A collection of base types for the DIFFEQ library.
Attempts an integration step for a single-step integrator.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(single_step_integrator), | intent(inout) | :: | this |
The single_step_integrator object. |
||
class(ode_container), | intent(inout) | :: | sys |
The ode_container object containing the ODE's to integrate. |
||
real(kind=real64), | intent(in) | :: | h |
The current step size. |
||
real(kind=real64), | intent(in) | :: | x |
The current value of the independent variable. |
||
real(kind=real64), | intent(in), | dimension(:) | :: | y |
An N-element array containing the current solution at x. |
|
real(kind=real64), | intent(in), | dimension(:) | :: | f |
An N-element array containing the values of the derivatives at x. |
|
real(kind=real64), | intent(out), | dimension(:) | :: | yn |
An N-element array where this routine will write the next solution estimate at x + h. |
|
real(kind=real64), | intent(out), | dimension(:) | :: | fn |
An N-element array where this routine will write the next derivative estimate at x + h. |
|
real(kind=real64), | intent(out), | dimension(:) | :: | yerr |
An N-element array where this routine will write an estimate of the error in each equation. |
|
real(kind=real64), | intent(out), | dimension(:,:) | :: | k |
An N-by-NSTAGES matrix containing the derivatives at each stage. |
Returns a logical parameter from a single_step_integrator object.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(single_step_integrator), | intent(in) | :: | this |
The single_step_integrator object. |
The parameter.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=real64), | intent(in) | :: | x | |||
real(kind=real64), | intent(in), | dimension(:) | :: | y | ||
real(kind=real64), | intent(out), | dimension(:) | :: | dydx |
Returns an integer value from the ode_integrator object.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(in) | :: | this |
The ode_integrator object. |
The requested value.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=real64), | intent(in) | :: | x | |||
real(kind=real64), | intent(in), | dimension(:) | :: | y | ||
real(kind=real64), | intent(out), | dimension(:,:) | :: | jac |
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=real64), | intent(in) | :: | x | |||
real(kind=real64), | intent(in), | dimension(:) | :: | y | ||
real(kind=real64), | intent(out), | dimension(:,:) | :: | m |
Solves the supplied system of ODE's.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ode_integrator), | intent(inout) | :: | this |
The ode_integrator object. |
||
class(ode_container), | intent(inout) | :: | sys |
The ode_container object containing the ODE's to integrate. |
||
real(kind=real64), | intent(in), | dimension(:) | :: | x |
An array, of at least 2 values, defining at a minimum the starting and ending values of the independent variable integration range. If more than two values are specified, the integration results will be returned at the supplied values. |
|
real(kind=real64), | intent(in), | dimension(:) | :: | iv |
An array containing the initial values for each ODE. |
|
class(errors), | intent(inout), | optional, | target | :: | err |
An optional errors-based object that if provided can be used to retrieve information relating to any errors encountered during execution. If not provided, a default implementation of the errors class is used internally to provide error handling. Possible errors and warning messages that may be encountered are as follows.
|
Gets an integer from the integrator.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(single_step_integrator), | intent(in) | :: | this |
The single_step_integrator object. |
The integer value.
Provides a routine for interpolation.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(single_step_integrator), | intent(in) | :: | this |
The single_step_integrator object. |
||
real(kind=real64), | intent(in) | :: | x |
The value of the independent variable at which to compute the interpolation. |
||
real(kind=real64), | intent(in) | :: | xn |
The previous value of the independent variable at which the solution is computed. |
||
real(kind=real64), | intent(in), | dimension(:) | :: | yn |
An N-element array containing the solution at xn. |
|
real(kind=real64), | intent(in), | dimension(:) | :: | fn |
An N-element array containing the derivatives at xn. |
|
real(kind=real64), | intent(in) | :: | xn1 |
The value of the independent variable at xn + h. |
||
real(kind=real64), | intent(in), | dimension(:) | :: | yn1 |
An N-element array containing the solution at xn + h. |
|
real(kind=real64), | intent(in), | dimension(:) | :: | fn1 |
An N-element array containing the derivatives at xn + h. |
|
real(kind=real64), | intent(out), | dimension(:) | :: | y |
An N-element array where this routine will write the solution values interpolated at x. |
Provides a routine for performing any actions, such as setting up interpolation, after successful completion of a step.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(single_step_integrator), | intent(inout) | :: | this |
The single_step_integrator object. |
||
class(ode_container), | intent(inout) | :: | sys |
The ode_container object containing the ODE's to integrate. |
||
logical, | intent(in) | :: | dense |
Determines if dense output is requested (true); else, false. |
||
real(kind=real64), | intent(in) | :: | x |
The previous value of the independent variable. |
||
real(kind=real64), | intent(in) | :: | xn |
The current value of the independent variable. |
||
real(kind=real64), | intent(in), | dimension(:) | :: | y |
An N-element array containing the solution at x. |
|
real(kind=real64), | intent(in), | dimension(:) | :: | yn |
An N-element array containing the solution at xn. |
|
real(kind=real64), | intent(in), | dimension(:) | :: | f |
An N-element array containing the derivatives at x. |
|
real(kind=real64), | intent(in), | dimension(:) | :: | fn |
An N-element array containing the derivatives at xn. |
|
real(kind=real64), | intent(inout), | dimension(:,:) | :: | k |
An N-by-NSTAGES matrix containing the derivatives at each stage. |
Provides a routine for performing any actions, such as setting up Jacobian calculations.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(single_step_integrator), | intent(inout) | :: | this |
The single_step_integrator object. |
||
logical, | intent(in) | :: | prevs |
Defines the status of the previous step. The value is true if the previous step was successful; else, false if the previous step failed. |
||
class(ode_container), | intent(inout) | :: | sys |
The ode_container object containing the ODE's to integrate. |
||
real(kind=real64), | intent(in) | :: | h |
The current step size. |
||
real(kind=real64), | intent(in) | :: | x |
The current value of the independent variable. |
||
real(kind=real64), | intent(in), | dimension(:) | :: | y |
An N-element array containing the current solution at x. |
|
real(kind=real64), | intent(in), | dimension(:) | :: | f |
An N-element array containing the values of the derivatives at x. |
|
class(errors), | intent(inout), | optional, | target | :: | err |
An optional errors-based object that if provided can be used to retrieve information relating to any errors encountered during execution. If not provided, a default implementation of the errors class is used internally to provide error handling. |
A container for the routine containing the ODEs to integrate.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
procedure(ode), | public, | pointer, nopass | :: | fcn | => | null() |
A pointer to the routine containing the ODEs to integrate. |
procedure(ode_jacobian), | public, | pointer, nopass | :: | jacobian | => | null() |
A pointer to the routine containing the analytical Jacobian. If supplied, this routine is utilized; however, if null, a finite difference approximation is utilized. |
procedure(ode_mass_matrix), | public, | pointer, nopass | :: | mass_matrix | => | null() |
A pointer to the routine containing the mass matrix for the system. If set to null (the default), an identity mass matrix will be assumed. |
procedure, public :: compute_jacobian => oc_jacobian | |
procedure, public :: get_finite_difference_step => oc_get_fd_step | |
procedure, public :: get_is_mass_matrix_dependent => oc_get_is_mass_dependent | |
procedure, public :: get_is_ode_defined => oc_get_is_ode_defined | |
procedure, public :: set_finite_difference_step => oc_set_fd_step | |
procedure, public :: set_is_mass_matrix_dependent => oc_set_is_mass_dependent |
The most basic ODE integrator object capable of integrating systems of ODE's.
procedure, public :: append_to_buffer => oi_append_to_buffer | ..\..\ Appends the supplied solution point to the internal solution buffer.<\p> |
procedure, public :: clear_buffer => oi_clear_buffer | ..\..\ Clears the contents of the buffer.<\p> |
procedure, public :: compute_error_norm => oi_estimate_error | ..\..\ Computes the norm of the scaled error estimate.<\p> |
procedure, public :: estimate_inital_step_size => oi_initial_step | ..\..\ Computes an estimate of an initial step size.<\p> |
procedure, public :: estimate_next_step_size => oi_next_step | ..\..\ Estimates the next step size.<\p> |
procedure, public :: get_absolute_tolerance => oi_get_abs_tol | ..\..\ Gets the absolute error tolerance.<\p> |
procedure, public :: get_allow_overshoot => oi_get_allow_overshoot | ..\..\ Gets a value determining if the solver is allowed to overshoot the final value in the integration range.<\p> |
procedure, public :: get_maximum_step_size => oi_get_max_step | ..\..\ Gets the magnitude of the maximum allowed step size.<\p> |
procedure, public :: get_minimum_step_size => oi_get_min_step | ..\..\ Gets the magnitude of the minimum allowed step size.<\p> |
procedure(ode_integer_inquiry), public, deferred, pass :: get_order | ..\..\ Returns the order of the integrator.<\p> |
procedure, public :: get_relative_tolerance => oi_get_abs_tol | ..\..\ Gets the relative error tolerance.<\p> |
procedure, public :: get_solution => oi_get_solution | ..\..\ Returns the solution computed by the integrator.<\p> |
procedure, public :: get_step_limit => oi_get_step_limit | ..\..\ Gets the limit on the number of integration steps.<\p> |
procedure, public :: get_step_size_control_parameter => oi_get_control_parameter | ..\..\ Gets the step size PI control parameter.<\p> |
procedure, public :: get_step_size_factor => oi_get_safety_factor | ..\..\ Gets the step size safety factor.<\p> |
procedure, public :: set_absolute_tolerance => oi_set_abs_tol | ..\..\ Sets the absolute error tolerance.<\p> |
procedure, public :: set_allow_overshoot => oi_set_allow_overshoot | ..\..\ Sets a value determining if the solver is allowed to overshoot the final value in the integration range.<\p> |
procedure, public :: set_maximum_step_size => oi_set_max_step | ..\..\ Sets the magnitude of the maximum allowed step size.<\p> |
procedure, public :: set_minimum_step_size => oi_set_min_step | ..\..\ Sets the magnitude of the minimum allowed step size.<\p> |
procedure, public :: set_relative_tolerance => oi_set_abs_tol | ..\..\ Sets the relative error tolerance.<\p> |
procedure, public :: set_step_limit => oi_set_step_limit | ..\..\ Sets the limit on the number of integration steps.<\p> |
procedure, public :: set_step_size_control_parameter => oi_set_control_parameter | ..\..\ Sets the step size PI control parameter.<\p> |
procedure, public :: set_step_size_factor => oi_set_safety_factor | ..\..\ Sets the step size safety factor.<\p> |
procedure(ode_solver), public, deferred, pass :: solve | ..\..\ Solves the supplied system of ODE's.<\p> |
The most basic, single-step integrator object capable of integrating systems of ODE's.
procedure, public :: append_to_buffer => oi_append_to_buffer | ..\..\ Appends the supplied solution point to the internal solution buffer.<\p> |
procedure(attempt_single_step), public, deferred, pass :: attempt_step | ..\..\ Attempts an integration step for a single-step integrator.<\p> |
procedure, public :: clear_buffer => oi_clear_buffer | ..\..\ Clears the contents of the buffer.<\p> |
procedure, public :: compute_error_norm => oi_estimate_error | ..\..\ Computes the norm of the scaled error estimate.<\p> |
procedure, public :: estimate_inital_step_size => oi_initial_step | ..\..\ Computes an estimate of an initial step size.<\p> |
procedure, public :: estimate_next_step_size => oi_next_step | ..\..\ Estimates the next step size.<\p> |
procedure, public :: get_absolute_tolerance => oi_get_abs_tol | ..\..\ Gets the absolute error tolerance.<\p> |
procedure, public :: get_allow_overshoot => oi_get_allow_overshoot | ..\..\ Gets a value determining if the solver is allowed to overshoot the final value in the integration range.<\p> |
procedure(get_single_step_logical_parameter), public, deferred, pass :: get_is_fsal | ..\..\ Gets a logical parameter stating if this is a first-same-as-last (FSAL) integrator.<\p> |
procedure, public :: get_maximum_step_size => oi_get_max_step | ..\..\ Gets the magnitude of the maximum allowed step size.<\p> |
procedure, public :: get_minimum_step_size => oi_get_min_step | ..\..\ Gets the magnitude of the minimum allowed step size.<\p> |
procedure(ode_integer_inquiry), public, deferred, pass :: get_order | ..\..\ Returns the order of the integrator.<\p> |
procedure, public :: get_relative_tolerance => oi_get_abs_tol | ..\..\ Gets the relative error tolerance.<\p> |
procedure, public :: get_solution => oi_get_solution | ..\..\ Returns the solution computed by the integrator.<\p> |
procedure(single_step_integer_inquiry), public, deferred, pass :: get_stage_count | ..\..\ Gets the number of stages used by the integrator.<\p> |
procedure, public :: get_step_limit => oi_get_step_limit | ..\..\ Gets the limit on the number of integration steps.<\p> |
procedure, public :: get_step_size_control_parameter => oi_get_control_parameter | ..\..\ Gets the step size PI control parameter.<\p> |
procedure, public :: get_step_size_factor => oi_get_safety_factor | ..\..\ Gets the step size safety factor.<\p> |
procedure(single_step_interpolate), public, deferred, pass :: interpolate | ..\..\ Performs an interpolation to estimate the solution at the requested point.<\p> |
procedure(single_step_post_step_routine), public, deferred, pass :: post_step_action | ..\..\ Performs actions such as setting up interpolation after completion of a successful integration step.<\p> |
procedure(single_step_pre_step_routine), public, deferred, pass :: pre_step_action | ..\..\ Provides a routine for performing any actions, such as setting up Jacobian calculations.<\p> |
procedure, public :: set_absolute_tolerance => oi_set_abs_tol | ..\..\ Sets the absolute error tolerance.<\p> |
procedure, public :: set_allow_overshoot => oi_set_allow_overshoot | ..\..\ Sets a value determining if the solver is allowed to overshoot the final value in the integration range.<\p> |
procedure, public :: set_maximum_step_size => oi_set_max_step | ..\..\ Sets the magnitude of the maximum allowed step size.<\p> |
procedure, public :: set_minimum_step_size => oi_set_min_step | ..\..\ Sets the magnitude of the minimum allowed step size.<\p> |
procedure, public :: set_relative_tolerance => oi_set_abs_tol | ..\..\ Sets the relative error tolerance.<\p> |
procedure, public :: set_step_limit => oi_set_step_limit | ..\..\ Sets the limit on the number of integration steps.<\p> |
procedure, public :: set_step_size_control_parameter => oi_set_control_parameter | ..\..\ Sets the step size PI control parameter.<\p> |
procedure, public :: set_step_size_factor => oi_set_safety_factor | ..\..\ Sets the step size safety factor.<\p> |
procedure, public :: solve => ssi_ode_solver | ..\..\ Solves the supplied system of ODE's.<\p> |